Contents
# # Contracting cube
#
# In this demo we simulate a unit cube that is fixed
# at $x = 0$ and free at $x = 1$. We use
# a transversally isotropic material with fiber oriented
# in the $x$-direction.
import dolfin
try:
    from dolfin_adjoint import (
        Constant,
        DirichletBC,
        Expression,
        UnitCubeMesh,
        interpolate,
        Mesh,
    )
except ImportError:
    from dolfin import (
        Constant,
        DirichletBC,
        interpolate,
        Expression,
        UnitCubeMesh,
        Mesh,
    )
import pulse
from fenics_plotly import plot
pulse.iterate.logger.setLevel(10)
# Create mesh
N = 6
mesh = UnitCubeMesh(N, N, N)
# Create subdomains
class Free(dolfin.SubDomain):
    def inside(self, x, on_boundary):
        return x[0] > (1.0 - dolfin.DOLFIN_EPS) and on_boundary


class Fixed(dolfin.SubDomain):
    def inside(self, x, on_boundary):
        return x[0] < dolfin.DOLFIN_EPS and on_boundary
# Create a facet fuction in order to mark the subdomains
ffun = dolfin.MeshFunction("size_t", mesh, 2)
ffun.set_all(0)
# Mark the first subdomain with value 1
fixed = Fixed()
fixed_marker = 1
fixed.mark(ffun, fixed_marker)
# Mark the second subdomain with value 2
free = Free()
free_marker = 2
free.mark(ffun, free_marker)
# Create a cell function (but we are not using it)
cfun = dolfin.MeshFunction("size_t", mesh, 3)
cfun.set_all(0)
# Collect the functions containing the markers
marker_functions = pulse.MarkerFunctions(ffun=ffun, cfun=cfun)
# Create mictrotructure
V_f = dolfin.VectorFunctionSpace(mesh, "CG", 1)
# Fibers
f0 = interpolate(Expression(("1.0", "0.0", "0.0"), degree=1), V_f)
# Sheets
s0 = interpolate(Expression(("0.0", "1.0", "0.0"), degree=1), V_f)
# Fiber-sheet normal
n0 = interpolate(Expression(("0.0", "0.0", "1.0"), degree=1), V_f)
# Collect the mictrotructure
microstructure = pulse.Microstructure(f0=f0, s0=s0, n0=n0)
# Create the geometry
geometry = pulse.Geometry(
    mesh=mesh,
    marker_functions=marker_functions,
    microstructure=microstructure,
)
# Use the default material parameters
material_parameters = pulse.HolzapfelOgden.default_parameters()
# Select model for active contraction
active_model = pulse.ActiveModels.active_strain
# active_model = "active_stress"
# Set the activation
activation = Constant(0.0)
# Create material
material = pulse.HolzapfelOgden(
    active_model=active_model,
    parameters=material_parameters,
    activation=activation,
)
# Make Dirichlet boundary conditions
def dirichlet_bc(W):
    V = W if W.sub(0).num_sub_spaces() == 0 else W.sub(0)
    return DirichletBC(V, Constant((0.0, 0.0, 0.0)), fixed)
# Make Neumann boundary conditions
neumann_bc = pulse.NeumannBC(traction=Constant(0.0), marker=free_marker)
# Collect Boundary Conditions
bcs = pulse.BoundaryConditions(dirichlet=(dirichlet_bc,), neumann=(neumann_bc,))
# Create problem
problem = pulse.MechanicsProblem(geometry, material, bcs)
# Solve problem
pulse.iterate.iterate(problem, activation, 0.1)
2021-11-23 07:02:17,244 [1335] DEBUG    pulse.iterate: Control: [0.0]
2021-11-23 07:02:17,245 [1335] DEBUG    pulse.iterate: Target: [0.1]
2021-11-23 07:02:17,247 [1335] DEBUG    pulse.iterate: Intial number of steps: 5 with step size 0.02
2021-11-23 07:02:17,248 [1335] INFO     pulse.iterate: Iterating to target control (f_21)...
2021-11-23 07:02:17,248 [1335] INFO     pulse.iterate: Current control: f_21 = 0.000
2021-11-23 07:02:17,249 [1335] INFO     pulse.iterate: Target: 0.100
2021-11-23 07:02:17,249 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:02:17,250 [1335] DEBUG    pulse.iterate: Maximum difference: 1.000e-01
2021-11-23 07:02:17,251 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:02:17,251 [1335] DEBUG    pulse.iterate: Current control: 0.020
2021-11-23 07:03:10,669 [1335] DEBUG    pulse.iterate: Solver did not converge...
2021-11-23 07:03:10,670 [1335] DEBUG    pulse.iterate: 
NOT CONVERGING
2021-11-23 07:03:10,671 [1335] DEBUG    pulse.iterate: Reduce control step
2021-11-23 07:03:10,672 [1335] DEBUG    pulse.iterate: Assign old state
2021-11-23 07:03:10,681 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:10,682 [1335] DEBUG    pulse.iterate: Maximum difference: 1.000e-01
2021-11-23 07:03:10,683 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:10,684 [1335] DEBUG    pulse.iterate: Current control: 0.010
2021-11-23 07:03:17,079 [1335] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2021-11-23 07:03:17,080 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:17,081 [1335] DEBUG    pulse.iterate: Maximum difference: 9.000e-02
2021-11-23 07:03:17,084 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:17,084 [1335] DEBUG    pulse.iterate: Maximum difference: 9.000e-02
2021-11-23 07:03:17,085 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:17,088 [1335] DEBUG    pulse.iterate: Current control: 0.020
2021-11-23 07:03:19,749 [1335] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2021-11-23 07:03:19,750 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:19,751 [1335] DEBUG    pulse.iterate: Maximum difference: 8.000e-02
2021-11-23 07:03:19,752 [1335] DEBUG    pulse.iterate: Adapt step size. New step size: 0.015
2021-11-23 07:03:19,753 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:19,754 [1335] DEBUG    pulse.iterate: Maximum difference: 8.000e-02
2021-11-23 07:03:19,755 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:19,756 [1335] DEBUG    pulse.iterate: Current control: 0.035
2021-11-23 07:03:23,539 [1335] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2021-11-23 07:03:23,540 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:23,541 [1335] DEBUG    pulse.iterate: Maximum difference: 6.500e-02
2021-11-23 07:03:23,542 [1335] DEBUG    pulse.iterate: Adapt step size. New step size: 0.022
2021-11-23 07:03:23,544 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:23,544 [1335] DEBUG    pulse.iterate: Maximum difference: 6.500e-02
2021-11-23 07:03:23,545 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:23,546 [1335] DEBUG    pulse.iterate: Current control: 0.058
2021-11-23 07:03:26,785 [1335] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2021-11-23 07:03:26,786 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:26,787 [1335] DEBUG    pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 07:03:26,788 [1335] DEBUG    pulse.iterate: Adapt step size. New step size: 0.034
2021-11-23 07:03:26,789 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:26,789 [1335] DEBUG    pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 07:03:26,790 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:26,791 [1335] DEBUG    pulse.iterate: Current control: 0.091
2021-11-23 07:03:28,395 [1335] DEBUG    pulse.iterate: Solver did not converge...
2021-11-23 07:03:28,396 [1335] DEBUG    pulse.iterate: 
NOT CONVERGING
2021-11-23 07:03:28,397 [1335] DEBUG    pulse.iterate: Reduce control step
2021-11-23 07:03:28,398 [1335] DEBUG    pulse.iterate: Assign old state
2021-11-23 07:03:28,406 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:28,406 [1335] DEBUG    pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 07:03:28,408 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:28,409 [1335] DEBUG    pulse.iterate: Current control: 0.074
2021-11-23 07:03:31,610 [1335] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2021-11-23 07:03:31,612 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:31,613 [1335] DEBUG    pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 07:03:31,614 [1335] DEBUG    pulse.iterate: Adapt step size. New step size: 0.025
2021-11-23 07:03:31,616 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:31,616 [1335] DEBUG    pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 07:03:31,617 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:31,618 [1335] DEBUG    pulse.iterate: Current control: 0.100
2021-11-23 07:03:33,261 [1335] DEBUG    pulse.iterate: Solver did not converge...
2021-11-23 07:03:33,262 [1335] DEBUG    pulse.iterate: 
NOT CONVERGING
2021-11-23 07:03:33,264 [1335] DEBUG    pulse.iterate: Reduce control step
2021-11-23 07:03:33,264 [1335] DEBUG    pulse.iterate: Assign old state
2021-11-23 07:03:33,271 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:33,272 [1335] DEBUG    pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 07:03:33,273 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:33,276 [1335] DEBUG    pulse.iterate: Current control: 0.087
2021-11-23 07:03:35,873 [1335] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2021-11-23 07:03:35,874 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:35,877 [1335] DEBUG    pulse.iterate: Maximum difference: 1.297e-02
2021-11-23 07:03:35,878 [1335] DEBUG    pulse.iterate: Adapt step size. New step size: 0.019
2021-11-23 07:03:35,880 [1335] DEBUG    pulse.iterate: Check target reached: NO
2021-11-23 07:03:35,882 [1335] DEBUG    pulse.iterate: Maximum difference: 1.297e-02
2021-11-23 07:03:35,882 [1335] DEBUG    pulse.iterate: Change step size for final iteration
2021-11-23 07:03:35,885 [1335] DEBUG    pulse.iterate: Try new control
2021-11-23 07:03:35,885 [1335] DEBUG    pulse.iterate: Current control: 0.100
2021-11-23 07:03:38,090 [1335] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2021-11-23 07:03:38,091 [1335] DEBUG    pulse.iterate: Check target reached: YES!
2021-11-23 07:03:38,094 [1335] DEBUG    pulse.iterate: Check target reached: YES!
([Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 47),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 195),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 246),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 313),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 372),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 466),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 552),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 595)],
 [Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 45),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 193),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 244),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 311),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 370),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 464),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 550),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 593)])
# Get displacement and hydrostatic pressure
u, p = problem.state.split(deepcopy=True)
V = dolfin.VectorFunctionSpace(mesh, "CG", 1)
u_int = interpolate(u, V)
new_mesh = Mesh(mesh)
dolfin.ALE.move(new_mesh, u_int)
fig = plot(mesh, show=False)
fig.add_plot(plot(new_mesh, color="red", show=False))
fig.show()